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

#include <G4NeutronHPFissionData.hh>

Inheritance diagram for G4NeutronHPFissionData:
G4VCrossSectionDataSet

Public Member Functions

 G4NeutronHPFissionData ()
 
 ~G4NeutronHPFissionData ()
 
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)
 
- 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 47 of file G4NeutronHPFissionData.hh.

Constructor & Destructor Documentation

G4NeutronHPFissionData::G4NeutronHPFissionData ( )

Definition at line 44 of file G4NeutronHPFissionData.cc.

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

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

Definition at line 59 of file G4NeutronHPFissionData.cc.

References G4PhysicsTable::clearAndDestroy().

60 {
61  if ( theCrossSections != NULL ) theCrossSections->clearAndDestroy();
62  delete theCrossSections;
63 }
void clearAndDestroy()

Member Function Documentation

void G4NeutronHPFissionData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 104 of file G4NeutronHPFissionData.cc.

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

105 {
106  if(&aP!=G4Neutron::Neutron())
107  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
108  size_t numberOfElements = G4Element::GetNumberOfElements();
109  //theCrossSections = new G4PhysicsTable( numberOfElements );
110  // TKDB
111  //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
112  if ( theCrossSections == NULL )
113  theCrossSections = new G4PhysicsTable( numberOfElements );
114  else
115  theCrossSections->clearAndDestroy();
116 
117  // make a PhysicsVector for each element
118 
119  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
120  for( size_t i=0; i<numberOfElements; ++i )
121  {
123  Instance()->MakePhysicsVector((*theElementTable)[i], this);
124  theCrossSections->push_back(physVec);
125  }
126 }
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:52
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NeutronHPData * Instance()
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
void clearAndDestroy()
void G4NeutronHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 128 of file G4NeutronHPFissionData.cc.

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

129 {
130  if(&aP!=G4Neutron::Neutron())
131  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
132 
133 //
134 // Dump element based cross section
135 // range 10e-5 eV to 20 MeV
136 // 10 point per decade
137 // in barn
138 //
139 
140  G4cout << G4endl;
141  G4cout << G4endl;
142  G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
143  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
144  G4cout << G4endl;
145  G4cout << "Name of Element" << G4endl;
146  G4cout << "Energy[eV] XS[barn]" << G4endl;
147  G4cout << G4endl;
148 
149  size_t numberOfElements = G4Element::GetNumberOfElements();
150  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
151 
152  for ( size_t i = 0 ; i < numberOfElements ; ++i )
153  {
154 
155  G4cout << (*theElementTable)[i]->GetName() << G4endl;
156 
157  if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
158  {
159  G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
160  G4cout << G4endl;
161  continue;
162  }
163 
164  G4int ie = 0;
165 
166  for ( ie = 0 ; ie < 130 ; ie++ )
167  {
168  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
169  G4bool outOfRange = false;
170 
171  if ( eKinetic < 20*MeV )
172  {
173  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
174  }
175 
176  }
177 
178  G4cout << G4endl;
179  }
180 
181  //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
182 }
#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
G4double G4NeutronHPFissionData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 187 of file G4NeutronHPFissionData.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::kelvin, G4ReactionProduct::Lorentz(), CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4Neutron::Neutron(), and G4ReactionProduct::SetMomentum().

Referenced by GetIsoCrossSection().

188 {
189  G4double result = 0;
190  if(anE->GetZ()<90) return result;
191  G4bool outOfRange;
192  G4int index = anE->GetIndex();
193 
194 // 100729 TK add safety
195 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
196 
197  // prepare neutron
198  G4double eKinetic = aP->GetKineticEnergy();
199  G4ReactionProduct theNeutron( aP->GetDefinition() );
200  theNeutron.SetMomentum( aP->GetMomentum() );
201  theNeutron.SetKineticEnergy( eKinetic );
202 
203  // prepare thermal nucleus
204  G4Nucleus aNuc;
205  G4double eps = 0.0001;
206  G4double theA = anE->GetN();
207  G4double theZ = anE->GetZ();
208  G4double eleMass;
209  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
211 
212  G4ReactionProduct boosted;
213  G4double aXsection;
214 
215  // MC integration loop
216  G4int counter = 0;
217  G4double buffer = 0;
218  G4int size = G4int(std::max(10., aT/60*kelvin));
219  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
220  G4double neutronVMag = neutronVelocity.mag();
221 
222  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
223  {
224  if(counter) buffer = result/counter;
225  while (counter<size)
226  {
227  counter ++;
228  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
229  boosted.Lorentz(theNeutron, aThermalNuc);
230  G4double theEkin = boosted.GetKineticEnergy();
231  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
232  // velocity correction.
233  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
234  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
235  result += aXsection;
236  }
237  size += size;
238  }
239  result /= counter;
240  return result;
241 }
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
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 G4NeutronHPFissionData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 78 of file G4NeutronHPFissionData.cc.

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

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 243 of file G4NeutronHPFissionData.cc.

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

244 {
246 }
static G4NeutronHPManager * GetInstance()
G4bool G4NeutronHPFissionData::IsIsoApplicable ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 65 of file G4NeutronHPFissionData.cc.

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

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 247 of file G4NeutronHPFissionData.cc.

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

248 {
250 }
static G4NeutronHPManager * GetInstance()
void SetVerboseLevel(G4int i)

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