Geant4  10.02.p03
G4ParticleHPFissionData Class Reference

#include <G4ParticleHPFissionData.hh>

Inheritance diagram for G4ParticleHPFissionData:
Collaboration diagram for G4ParticleHPFissionData:

Public Member Functions

 G4ParticleHPFissionData ()
 
 ~G4ParticleHPFissionData ()
 
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 &)
 
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 48 of file G4ParticleHPFissionData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFissionData()

G4ParticleHPFissionData::G4ParticleHPFissionData ( )

Definition at line 47 of file G4ParticleHPFissionData.cc.

48 :G4VCrossSectionDataSet("NeutronHPFissionXS")
49 {
50  SetMinKinEnergy( 0*MeV );
51  SetMaxKinEnergy( 20*MeV );
52 
53  theCrossSections = 0;
54  onFlightDB = true;
55  instanceOfWorker = false;
57  instanceOfWorker = true;
58  }
59  //BuildPhysicsTable(*G4Neutron::Neutron());
60 }
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:

◆ ~G4ParticleHPFissionData()

G4ParticleHPFissionData::~G4ParticleHPFissionData ( )

Definition at line 62 of file G4ParticleHPFissionData.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPFissionData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 110 of file G4ParticleHPFissionData.cc.

111 {
112 
113  if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() ) {
114  G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
115  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of fission reaction of neutrons (<20MeV)." << G4endl;
116  onFlightDB = false;
117  }
118 
119  if(&aP!=G4Neutron::Neutron())
120  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
121 
122  if ( G4Threading::IsWorkerThread() ) {
124  return;
125  }
126 
127  size_t numberOfElements = G4Element::GetNumberOfElements();
128  //theCrossSections = new G4PhysicsTable( numberOfElements );
129  // TKDB
130  //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
131  if ( theCrossSections == NULL )
132  theCrossSections = new G4PhysicsTable( numberOfElements );
133  else
135 
136  // make a PhysicsVector for each element
137 
138  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
139  for( size_t i=0; i<numberOfElements; ++i )
140  {
142  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
143  theCrossSections->push_back(physVec);
144  }
145 
147 }
static G4ParticleHPManager * GetInstance()
void RegisterFissionCrossSections(G4PhysicsTable *val)
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:89
G4PhysicsTable * GetFissionCrossSections()
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
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 G4ParticleHPFissionData::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 284 of file G4ParticleHPFissionData.cc.

285 {
286  outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n" ;
287 }

◆ DumpPhysicsTable()

void G4ParticleHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 149 of file G4ParticleHPFissionData.cc.

150 {
151  if(&aP!=G4Neutron::Neutron())
152  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
153 
154 //
155 // Dump element based cross section
156 // range 10e-5 eV to 20 MeV
157 // 10 point per decade
158 // in barn
159 //
160 
161  G4cout << G4endl;
162  G4cout << G4endl;
163  G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
164  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
165  G4cout << G4endl;
166  G4cout << "Name of Element" << G4endl;
167  G4cout << "Energy[eV] XS[barn]" << G4endl;
168  G4cout << G4endl;
169 
170  size_t numberOfElements = G4Element::GetNumberOfElements();
171  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
172 
173  for ( size_t i = 0 ; i < numberOfElements ; ++i )
174  {
175 
176  G4cout << (*theElementTable)[i]->GetName() << G4endl;
177 
178  if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
179  {
180  G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
181  G4cout << G4endl;
182  continue;
183  }
184 
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  }
198 
199  G4cout << G4endl;
200  }
201 
202  //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
203 }
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:

◆ GetCrossSection()

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

Definition at line 208 of file G4ParticleHPFissionData.cc.

209 {
210  G4double result = 0;
211  if(anE->GetZ()<88) return result;
212  G4bool outOfRange;
213  G4int index = anE->GetIndex();
214 
215 // 100729 TK add safety
216 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
217 
218  // prepare neutron
219  G4double eKinetic = aP->GetKineticEnergy();
220  G4ReactionProduct theNeutronRP( aP->GetDefinition() );
221  theNeutronRP.SetMomentum( aP->GetMomentum() );
222  theNeutronRP.SetKineticEnergy( eKinetic );
223 
224  if ( !onFlightDB ) {
225  //NEGLECT_DOPPLER
226  G4double factor = 1.0;
227  if ( eKinetic < aT * k_Boltzmann ) {
228  // below 0.1 eV neutrons
229  // Have to do some, but now just igonre.
230  // Will take care after performance check.
231  // factor = factor * targetV;
232  }
233  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
234  }
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  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
244 
245  G4ReactionProduct boosted;
246  G4double aXsection;
247 
248  // MC integration loop
249  G4int counter = 0;
250  G4double buffer = 0;
251  G4int size = G4int(std::max(10., aT/60*kelvin));
252  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
253  G4double neutronVMag = neutronVelocity.mag();
254 
255  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
256  {
257  if(counter) buffer = result/counter;
258  while (counter<size) // Loop checking, 11.05.2015, T. Koi
259  {
260  counter ++;
261  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
262  boosted.Lorentz(theNeutronRP, aThermalNuc);
263  G4double theEkin = boosted.GetKineticEnergy();
264  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
265  // velocity correction.
266  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
267  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
268  result += aXsection;
269  }
270  size += size;
271  }
272  result /= counter;
273  return result;
274 }
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 G4ParticleHPFissionData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 84 of file G4ParticleHPFissionData.cc.

89 {
90  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
91 
92  ke_cache = dp->GetKineticEnergy();
93  element_cache = element;
95  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
96  xs_cache = xs;
97  return xs;
98 }
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 G4ParticleHPFissionData::GetVerboseLevel ( ) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 276 of file G4ParticleHPFissionData.cc.

277 {
279 }
static G4ParticleHPManager * GetInstance()
Here is the call graph for this function:

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 71 of file G4ParticleHPFissionData.cc.

75 {
76  G4double eKin = dp->GetKineticEnergy();
77  if ( eKin > GetMaxKinEnergy()
78  || eKin < GetMinKinEnergy()
79  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
80 
81  return true;
82 }
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 G4ParticleHPFissionData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 280 of file G4ParticleHPFissionData.cc.

281 {
283 }
static G4ParticleHPManager * GetInstance()
Here is the call graph for this function:

Member Data Documentation

◆ element_cache

const G4Element* G4ParticleHPFissionData::element_cache
private

Definition at line 93 of file G4ParticleHPFissionData.hh.

◆ instanceOfWorker

G4bool G4ParticleHPFissionData::instanceOfWorker
private

Definition at line 89 of file G4ParticleHPFissionData.hh.

◆ ke_cache

G4double G4ParticleHPFissionData::ke_cache
private

Definition at line 91 of file G4ParticleHPFissionData.hh.

◆ material_cache

const G4Material* G4ParticleHPFissionData::material_cache
private

Definition at line 94 of file G4ParticleHPFissionData.hh.

◆ onFlightDB

G4bool G4ParticleHPFissionData::onFlightDB
private

Definition at line 88 of file G4ParticleHPFissionData.hh.

◆ theCrossSections

G4PhysicsTable* G4ParticleHPFissionData::theCrossSections
private

Definition at line 86 of file G4ParticleHPFissionData.hh.

◆ xs_cache

G4double G4ParticleHPFissionData::xs_cache
private

Definition at line 92 of file G4ParticleHPFissionData.hh.


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